#================================================
# A km-sub model constructed using MultilayerPy (pre public release)
# MultilayerPy v1.0 can be accessed here (https://doi.org/10.5281/zenodo.6411189)
# A separate Jupyter Notebook which will run the model code below will be released with v1.1 in the source code (this will be available via the above doi)
#
# **How to run manually using scipy.integrate.solve_ivp:**
#
# >>> from scipy.integrate import solve_ivp
#
# param_dict: a dictionary with relevant parameters as keys and their values as values (see below)
# V: array with the volume of each model bulk layer
# A: array with surface area of each model bulk layer
# Lorg: The number of model bulk layers
# layer_thick: array of the thickness of each model bulk layer
# 
# y0:array of initial surface/bulk concentrations of each component of length Lorg * 7 + 7 (there are 7 components) in component order: 
#    surface concs. (cm-2), bulk concs. (cm-3). 
#    oleic acid initial surf. conc. = 9.68e13 cm-2, oleic acid initial bulk conc. = 1.21e21 cm-3, the rest are 0
#
# >>> y0 = [oleic acid initital surf conc, Lorg * oleic acid initial bulk conc., ***Repeat for other 6 components***]
#
# t_span: a tuple with t0 (start time of integration) and tf (end time of integration), in seconds
#
# >>> t_span = (t0,tf)
#
# >>> args = (param_dict,V,A,Lorg,layer_thick)
# >>> result = solve_ivp(dydt,t_span,y0,method='BDF',args=args)	
#	
# consult the scipy.integrate.solve_ivp documentation (https://docs.scipy.org/doc/scipy/reference/generated/scipy.integrate.solve_ivp.html)
# for full instructions on how to access the model output. Summary:
#
# result.t --> array of saved time points for model output
# result.y --> array of shape (Lorg * 7 + 7, len(result.t)), surface and bulk concentrations in each layer at each time point
#
# Reaction name: Oleic acid ozonolysis lamellar
# Geometry: film
# Number of model components: 7
# Diffusion regime: vignes
#================================================

import numpy as np
def dydt(t,y,param_dict,V,A,Lorg,layer_thick,param_evolution_func=None,additional_params=None):
    """ Function defining ODEs, returns list of dydt"""

    if typec!= type(None):
        param_dict = param_evolution_func(t,y,param_dict,additional_params=additional_params)

    #--------------Unpack parameters (random order)---------------

    # If MultilayerPy used, this "try" block works, otherwise it will run the "except" block
    # see except block for parameter descriptions
    try: 
        delta_3 = param_dict["delta_3"].value
        delta_4 = param_dict["delta_4"].value
        k_1_2 = param_dict["k_1_2"].value
        k_1_2_surf = param_dict["k_1_2_surf"].value
        Db_5 = param_dict["Db_5"].value
        Db_5_7 = param_dict["Db_5_7"].value
        Db_7 = param_dict["Db_7"].value
        Db_3 = param_dict["Db_3"].value
        w_2 = param_dict["w_2"].value
        delta_5 = param_dict["delta_5"].value
        Db_1_7 = param_dict["Db_1_7"].value
        Db_4_7 = param_dict["Db_4_7"].value
        Db_3_7 = param_dict["Db_3_7"].value
        Db_2 = param_dict["Db_2"].value
        Db_2_7 = param_dict["Db_2_7"].value
        delta_2 = param_dict["delta_2"].value
        Db_1_6 = param_dict["Db_1_6"].value
        Td_2 = param_dict["Td_2"].value
        delta_6 = param_dict["delta_6"].value
        Db_3_6 = param_dict["Db_3_6"].value
        Db_6 = param_dict["Db_6"].value
        Xgs_2 = param_dict["Xgs_2"].value
        k_3_6 = param_dict["k_3_6"].value
        k_3_6_surf = param_dict["k_3_6_surf"].value
        k_3_4 = param_dict["k_3_4"].value
        k_3_4_surf = param_dict["k_3_4_surf"].value
        Db_1 = param_dict["Db_1"].value
        Db_2_6 = param_dict["Db_2_6"].value
        T = param_dict["T"].value
        Db_4 = param_dict["Db_4"].value
        H_2 = param_dict["H_2"].value
        delta_7 = param_dict["delta_7"].value
        alpha_s_0_2 = param_dict["alpha_s_0_2"].value
        Db_4_6 = param_dict["Db_4_6"].value
        delta_1 = param_dict["delta_1"].value
        Db_5_6 = param_dict["Db_5_6"].value

    except AttributeError:
        delta_3 = param_dict["delta_3"] # molec. diameter of CI (cm)
        delta_4 = param_dict["delta_4"] # molec. diameter of C9 (cm)
        k_1_2 = param_dict["k_1_2"] # bulk rxn constant - oleic acid + o3 (cm3 s-1)
        k_1_2_surf = param_dict["k_1_2_surf"] # surf rxn constant - oleic acid + o3 (cm2 s-1)
        Db_5 = param_dict["Db_5"] # bulk diffusion coeff. nonanal (cm2 s-1)
        Db_5_7 = param_dict["Db_5_7"] # bulk diffusion coeff nonanal in trimer (cm2 s-1)
        Db_7 = param_dict["Db_7"] # bulk diffusion coeff trimer (cm2 s-1)
        Db_3 = param_dict["Db_3"] # bulk diffusion coeff CI (cm2 s-1)
        w_2 = param_dict["w_2"] # mean thermal velocity O3 (cm s-1)
        delta_5 = param_dict["delta_5"] # molec. diameter of nonanal (cm)
        Db_1_7 = param_dict["Db_1_7"] # bulk diffusion coeff of oleic acid in trimer (cm2 s-1)
        Db_4_7 = param_dict["Db_4_7"] # bulk diffusion coeff of C9 in trimer (cm2 s-1)
        Db_3_7 = param_dict["Db_3_7"] # bulk diffusion coeff of CI in trimer (cm2 s-1)
        Db_2 = param_dict["Db_2"] # bulk diffusion coeff of O3
        Db_2_7 = param_dict["Db_2_7"] # bulk diffusion coeff of O3 in trimer (cm2 s-1)
        delta_2 = param_dict["delta_2"] # molec. diameter of O3 (cm)
        Db_1_6 = param_dict["Db_1_6"] # bulk diffusion coeff of oleic acid in dimer (cm2 s-1)
        Td_2 = param_dict["Td_2"] # surface desorption lifetime of O3
        delta_6 = param_dict["delta_6"] # molec. diameter of dimer
        Db_3_6 = param_dict["Db_3_6"] # bulk diffusion coeff of CI in dimer (cm2 s-1)
        Db_6 = param_dict["Db_6"] # bulk diffusion coeff of dimer (cm2 s-1)
        Xgs_2 = param_dict["Xgs_2"] # near surface gas phase conc. of O3 (cm-3)
        k_3_6 = param_dict["k_3_6"] # bulk rxn constant - CI + dimer (cm3 s-1)
        k_3_6_surf = param_dict["k_3_6_surf"] # surf rxn constant - CI + dimer (cm2 s-1)
        k_3_4 = param_dict["k_3_4"] # bulk rxn constant - CI + C9 (cm3 s-1)
        k_3_4_surf = param_dict["k_3_4_surf"] # surf rxn constant - CI + C9 (cm2 s-1)
        Db_1 = param_dict["Db_1"] # bulk diffusion coeff oleic acid (cm2 s-1)
        Db_2_6 = param_dict["Db_2_6"] # bulk diffusion coeff O3 in dimer (cm2 s-1)
        T = param_dict["T"] # Temp. (K)
        Db_4 = param_dict["Db_4"] # bulk diffusion coeff C9 (cm2 s-1)
        H_2 = param_dict["H_2"] # Henry's law coeff O3 (mol cm-3 atm-1)
        delta_7 = param_dict["delta_7"] # molec. diameter trimer (cm)
        alpha_s_0_2 = param_dict["alpha_s_0_2"] # surface accom. coeff O3 on clear surface
        Db_4_6 = param_dict["Db_4_6"] # bulk diffusion coeff C9 in dimer (cm2 s-1)
        delta_1 = param_dict["delta_1"] # molec. dimeter oleic acid (cm)
        Db_5_6 = param_dict["Db_5_6"] # bulk diffusion coeff nonanal in dimer (cm2 s-1)

    R = 82.0578
    Na = 6.022e23

    # init dydt array
    dydt = np.zeros(Lorg * 7 + 7)

    #--------------Define surface uptake parameters for gas components---------------

    #calculate surface fraction of each component
    fs_1 = y[0*Lorg+0] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_2 = y[1*Lorg+1] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_3 = y[2*Lorg+2] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_4 = y[3*Lorg+3] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_5 = y[4*Lorg+4] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_6 = y[5*Lorg+5] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])
    fs_7 = y[6*Lorg+6] / (y[0*Lorg+0] + y[1*Lorg+1]+ y[2*Lorg+2]+ y[3*Lorg+3]+ y[4*Lorg+4]+ y[5*Lorg+5]+ y[6*Lorg+6])


    # component 2 surf params

    surf_cover = delta_2**2 * y[1*Lorg+1] 
    alpha_s_2 = alpha_s_0_2 * (1-surf_cover)
    J_coll_X_2 = Xgs_2 * w_2/4
    J_ads_X_2 = alpha_s_2 * J_coll_X_2
    J_des_X_2 = (1 / Td_2) * y[1*Lorg+1]
    #--------------Bulk Diffusion evolution---------------

    # Db and fb arrays
    fb_1 = y[0*Lorg+1:1*Lorg+0+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_2 = y[1*Lorg+2:2*Lorg+1+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_3 = y[2*Lorg+3:3*Lorg+2+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_4 = y[3*Lorg+4:4*Lorg+3+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_5 = y[4*Lorg+5:5*Lorg+4+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_6 = y[5*Lorg+6:6*Lorg+5+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])
    fb_7 = y[6*Lorg+7:7*Lorg+6+1] / (y[0*Lorg+1:1*Lorg+0+1] + y[1*Lorg+2:2*Lorg+1+1]+ y[2*Lorg+3:3*Lorg+2+1]+ y[3*Lorg+4:4*Lorg+3+1]+ y[4*Lorg+5:5*Lorg+4+1]+ y[5*Lorg+6:6*Lorg+5+1]+ y[6*Lorg+7:7*Lorg+6+1])

    Db_1_arr = np.ones(Lorg) * Db_1
    Db_2_arr = np.ones(Lorg) * Db_2
    Db_3_arr = np.ones(Lorg) * Db_3
    Db_4_arr = np.ones(Lorg) * Db_4
    Db_5_arr = np.ones(Lorg) * Db_5
    Db_6_arr = np.ones(Lorg) * Db_6
    Db_7_arr = np.ones(Lorg) * Db_7
    Db_5_7_arr = np.ones(Lorg) * Db_5_7
    Db_1_7_arr = np.ones(Lorg) * Db_1_7
    Db_4_7_arr = np.ones(Lorg) * Db_4_7
    Db_3_7_arr = np.ones(Lorg) * Db_3_7
    Db_2_7_arr = np.ones(Lorg) * Db_2_7
    Db_1_6_arr = np.ones(Lorg) * Db_1_6
    Db_3_6_arr = np.ones(Lorg) * Db_3_6
    Db_2_6_arr = np.ones(Lorg) * Db_2_6
    Db_4_6_arr = np.ones(Lorg) * Db_4_6
    Db_5_6_arr = np.ones(Lorg) * Db_5_6

    # bulk diffusion

    Db_1_arr = (Db_1_arr**(1-fb_6-fb_7)) * (Db_1_6_arr**fb_6) * (Db_1_7_arr**fb_7) 
    Db_2_arr = (Db_2_arr**(1-fb_6-fb_7)) * (Db_2_6_arr**fb_6) * (Db_2_7_arr**fb_7) 
    Db_3_arr = (Db_3_arr**(1-fb_6-fb_7)) * (Db_3_6_arr**fb_6) * (Db_3_7_arr**fb_7) 
    Db_4_arr = (Db_4_arr**(1-fb_6-fb_7)) * (Db_4_6_arr**fb_6) * (Db_4_7_arr**fb_7) 
    Db_5_arr = (Db_5_arr**(1-fb_6-fb_7)) * (Db_5_6_arr**fb_6) * (Db_5_7_arr**fb_7) 
    Db_6_arr = Db_6_arr
    Db_7_arr = Db_7_arr
    kbby_1 = (4/np.pi) * Db_1_arr / layer_thick 
    kbbx_2 = (4/np.pi) * Db_2_arr / layer_thick 
    kbby_3 = (4/np.pi) * Db_3_arr / layer_thick 
    kbby_4 = (4/np.pi) * Db_4_arr / layer_thick 
    kbby_5 = (4/np.pi) * Db_5_arr / layer_thick 
    kbby_6 = (4/np.pi) * Db_6_arr / layer_thick 
    kbby_7 = (4/np.pi) * Db_7_arr / layer_thick 


    # surface diffusion
    Ds_1 = (Db_1**fs_1) * (Db_1_6**fs_6) * (Db_1_7**fs_7) 
    Ds_2 = (Db_2**fs_2) * (Db_2_6**fs_6) * (Db_2_7**fs_7) 
    Ds_3 = (Db_3**fs_3) * (Db_3_6**fs_6) * (Db_3_7**fs_7) 
    Ds_4 = (Db_4**fs_4) * (Db_4_6**fs_6) * (Db_4_7**fs_7) 
    Ds_5 = (Db_5**fs_5) * (Db_5_6**fs_6) * (Db_5_7**fs_7) 
    Ds_6 = Db_6
    Ds_7 = Db_7
    kbs_2 = (8/np.pi) * Ds_2 / (layer_thick[0] + delta_1+ delta_2+ delta_3+ delta_4+ delta_5+ delta_6+ delta_7)
    H_cc_2 = H_2 * R * T
    ksb_2 = H_cc_2 * kbs_2 / Td_2 / (w_2*alpha_s_2/4) 
    kbss_1 = (8*Db_1_arr[0])/((layer_thick[0]+delta_1)*np.pi) 
    kbss_3 = (8*Db_3_arr[0])/((layer_thick[0]+delta_3)*np.pi) 
    kbss_4 = (8*Db_4_arr[0])/((layer_thick[0]+delta_4)*np.pi) 
    kbss_5 = (8*Db_5_arr[0])/((layer_thick[0]+delta_5)*np.pi) 
    kbss_6 = (8*Db_6_arr[0])/((layer_thick[0]+delta_6)*np.pi) 
    kbss_7 = (8*Db_7_arr[0])/((layer_thick[0]+delta_7)*np.pi) 
    kssb_1 = kbss_1 / delta_1 
    kssb_3 = kbss_3 / delta_3 
    kssb_4 = kbss_4 / delta_4 
    kssb_5 = kbss_5 / delta_5 
    kssb_6 = kbss_6 / delta_6 
    kssb_7 = kbss_7 / delta_7 


    # mass fluxes

    Jssb_1 = kssb_1 * y[0*Lorg+0] 
    Jbss_1 = kbss_1 * y[0*Lorg+1] 
    Jb1b2_1 = kbby_1[0] * (y[0*Lorg+1+1] - y[0*Lorg+1]) 
    Jbb_core_1 = kbby_1[-1] * (y[1*Lorg+0-1] - y[1*Lorg+0]) 
    Jsb_2 = ksb_2 * y[1*Lorg+1] 
    Jbs_2 = kbs_2 * y[1*Lorg+2] 
    Jb1b2_2 = kbbx_2[0] * (y[1*Lorg+2+1] - y[1*Lorg+2]) 
    Jbb_core_2 = kbbx_2[-1] * (y[2*Lorg+1-1] - y[2*Lorg+1]) 
    Jssb_3 = kssb_3 * y[2*Lorg+2] 
    Jbss_3 = kbss_3 * y[2*Lorg+3] 
    Jb1b2_3 = kbby_3[0] * (y[2*Lorg+3+1] - y[2*Lorg+3]) 
    Jbb_core_3 = kbby_3[-1] * (y[3*Lorg+2-1] - y[3*Lorg+2]) 
    Jssb_4 = kssb_4 * y[3*Lorg+3] 
    Jbss_4 = kbss_4 * y[3*Lorg+4] 
    Jb1b2_4 = kbby_4[0] * (y[3*Lorg+4+1] - y[3*Lorg+4]) 
    Jbb_core_4 = kbby_4[-1] * (y[4*Lorg+3-1] - y[4*Lorg+3]) 
    Jssb_5 = kssb_5 * y[4*Lorg+4] 
    Jbss_5 = kbss_5 * y[4*Lorg+5] 
    Jb1b2_5 = kbby_5[0] * (y[4*Lorg+5+1] - y[4*Lorg+5]) 
    Jbb_core_5 = kbby_5[-1] * (y[5*Lorg+4-1] - y[5*Lorg+4]) 
    Jssb_6 = kssb_6 * y[5*Lorg+5] 
    Jbss_6 = kbss_6 * y[5*Lorg+6] 
    Jb1b2_6 = kbby_6[0] * (y[5*Lorg+6+1] - y[5*Lorg+6]) 
    Jbb_core_6 = kbby_6[-1] * (y[6*Lorg+5-1] - y[6*Lorg+5]) 
    Jssb_7 = kssb_7 * y[6*Lorg+6] 
    Jbss_7 = kbss_7 * y[6*Lorg+7] 
    Jb1b2_7 = kbby_7[0] * (y[6*Lorg+7+1] - y[6*Lorg+7]) 
    Jbb_core_7 = kbby_7[-1] * (y[7*Lorg+6-1] - y[7*Lorg+6]) 

    c = 0.454

    #----component number 1, Oleic acid----
    dydt[0*Lorg+0] = Jbss_1 - Jssb_1 - y[0*Lorg+0] * y[1*Lorg+1] * k_1_2_surf
    dydt[0*Lorg+1] = (Jssb_1 - Jbss_1) * (A[0]/V[0]) + Jb1b2_1 * (A[1]/V[0]) - y[0*Lorg+1] * y[1*Lorg+2] * k_1_2
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_1 = kbby_1[i] * (y[0*Lorg+1+i-1] - y[0*Lorg+1+i]) 
        Jbb_plus1_1 = kbby_1[i+1] * (y[0*Lorg+1+i+1] - y[0*Lorg+1+i]) 
        dydt[0*Lorg+1+i] = Jbb_minus1_1 * (A[i]/V[i]) + Jbb_plus1_1 * (A[i+1]/V[i]) - y[0*Lorg+1+i] * y[1*Lorg+2+i] * k_1_2
    dydt[1*Lorg+0] = Jbb_core_1 * (A[-1]/V[-1]) - y[1*Lorg+0] * y[2*Lorg+1] * k_1_2

    #----component number 2, Ozone----
    dydt[1*Lorg+1] = Jbs_2 - Jsb_2 + J_ads_X_2 - J_des_X_2- y[1*Lorg+1] * y[0*Lorg+0] * k_1_2_surf
    dydt[1*Lorg+2] = (Jsb_2 - Jbs_2) * (A[0]/V[0]) + Jb1b2_2 * (A[1]/V[0]) - y[1*Lorg+2] * y[0*Lorg+1] * k_1_2
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_2 = kbbx_2[i] * (y[1*Lorg+2+i-1] - y[1*Lorg+2+i]) 
        Jbb_plus1_2 = kbbx_2[i+1] * (y[1*Lorg+2+i+1] - y[1*Lorg+2+i]) 
        dydt[1*Lorg+2+i] = Jbb_minus1_2 * (A[i]/V[i]) + Jbb_plus1_2 * (A[i+1]/V[i]) - y[1*Lorg+2+i] * y[0*Lorg+1+i] * k_1_2
    dydt[2*Lorg+1] = Jbb_core_2 * (A[-1]/V[-1]) - y[2*Lorg+1] * y[1*Lorg+0] * k_1_2

    #----component number 3, CI----
    dydt[2*Lorg+2] = Jbss_3 - Jssb_3 + 1.0 * y[0*Lorg+0] * y[1*Lorg+1] * k_1_2_surf- y[2*Lorg+2] * y[3*Lorg+3] * k_3_4_surf- y[2*Lorg+2] * y[5*Lorg+5] * k_3_6_surf
    dydt[2*Lorg+3] = (Jssb_3 - Jbss_3) * (A[0]/V[0]) + Jb1b2_3 * (A[1]/V[0]) + 1.0 * y[0*Lorg+1] * y[1*Lorg+2] * k_1_2- y[2*Lorg+3] * y[3*Lorg+4] * k_3_4- y[2*Lorg+3] * y[5*Lorg+6] * k_3_6
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_3 = kbby_3[i] * (y[2*Lorg+3+i-1] - y[2*Lorg+3+i]) 
        Jbb_plus1_3 = kbby_3[i+1] * (y[2*Lorg+3+i+1] - y[2*Lorg+3+i]) 
        dydt[2*Lorg+3+i] = Jbb_minus1_3 * (A[i]/V[i]) + Jbb_plus1_3 * (A[i+1]/V[i]) + 1.0 * y[0*Lorg+1+i] * y[1*Lorg+2+i] * k_1_2- y[2*Lorg+3+i] * y[3*Lorg+4+i] * k_3_4- y[2*Lorg+3+i] * y[5*Lorg+6+i] * k_3_6
    dydt[3*Lorg+2] = Jbb_core_3 * (A[-1]/V[-1]) + 1.0 * y[1*Lorg+0] * y[2*Lorg+1] * k_1_2- y[3*Lorg+2] * y[4*Lorg+3] * k_3_4- y[3*Lorg+2] * y[6*Lorg+5] * k_3_6

    #----component number 4, C9----
    dydt[3*Lorg+3] = Jbss_4 - Jssb_4 + (1-c) * y[0*Lorg+0] * y[1*Lorg+1] * k_1_2_surf- y[3*Lorg+3] * y[2*Lorg+2] * k_3_4_surf
    dydt[3*Lorg+4] = (Jssb_4 - Jbss_4) * (A[0]/V[0]) + Jb1b2_4 * (A[1]/V[0]) + (1-c) * y[0*Lorg+1] * y[1*Lorg+2] * k_1_2- y[3*Lorg+4] * y[2*Lorg+3] * k_3_4
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_4 = kbby_4[i] * (y[3*Lorg+4+i-1] - y[3*Lorg+4+i]) 
        Jbb_plus1_4 = kbby_4[i+1] * (y[3*Lorg+4+i+1] - y[3*Lorg+4+i]) 
        dydt[3*Lorg+4+i] = Jbb_minus1_4 * (A[i]/V[i]) + Jbb_plus1_4 * (A[i+1]/V[i]) + (1-c) * y[0*Lorg+1+i] * y[1*Lorg+2+i] * k_1_2- y[3*Lorg+4+i] * y[2*Lorg+3+i] * k_3_4
    dydt[4*Lorg+3] = Jbb_core_4 * (A[-1]/V[-1]) + (1-c) * y[1*Lorg+0] * y[2*Lorg+1] * k_1_2- y[4*Lorg+3] * y[3*Lorg+2] * k_3_4

    kloss = 1.4
    #----component number 5, Nonanal----
    dydt[4*Lorg+4] = Jbss_5 - Jssb_5 + c * y[0*Lorg+0] * y[1*Lorg+1] * k_1_2_surf - kloss * y[4*Lorg+4]
    dydt[4*Lorg+5] = (Jssb_5 - Jbss_5) * (A[0]/V[0]) + Jb1b2_5 * (A[1]/V[0]) + c * y[0*Lorg+1] * y[1*Lorg+2] * k_1_2
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_5 = kbby_5[i] * (y[4*Lorg+5+i-1] - y[4*Lorg+5+i]) 
        Jbb_plus1_5 = kbby_5[i+1] * (y[4*Lorg+5+i+1] - y[4*Lorg+5+i]) 
        dydt[4*Lorg+5+i] = Jbb_minus1_5 * (A[i]/V[i]) + Jbb_plus1_5 * (A[i+1]/V[i]) + c * y[0*Lorg+1+i] * y[1*Lorg+2+i] * k_1_2
    dydt[5*Lorg+4] = Jbb_core_5 * (A[-1]/V[-1]) + c * y[1*Lorg+0] * y[2*Lorg+1] * k_1_2

    #----component number 6, Dimer----
    dydt[5*Lorg+5] = Jbss_6 - Jssb_6 + 1.0 * y[2*Lorg+2] * y[3*Lorg+3] * k_3_4_surf- y[5*Lorg+5] * y[2*Lorg+2] * k_3_6_surf
    dydt[5*Lorg+6] = (Jssb_6 - Jbss_6) * (A[0]/V[0]) + Jb1b2_6 * (A[1]/V[0]) + 1.0 * y[2*Lorg+3] * y[3*Lorg+4] * k_3_4- y[5*Lorg+6] * y[2*Lorg+3] * k_3_6
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_6 = kbby_6[i] * (y[5*Lorg+6+i-1] - y[5*Lorg+6+i]) 
        Jbb_plus1_6 = kbby_6[i+1] * (y[5*Lorg+6+i+1] - y[5*Lorg+6+i]) 
        dydt[5*Lorg+6+i] = Jbb_minus1_6 * (A[i]/V[i]) + Jbb_plus1_6 * (A[i+1]/V[i]) + 1.0 * y[2*Lorg+3+i] * y[3*Lorg+4+i] * k_3_4- y[5*Lorg+6+i] * y[2*Lorg+3+i] * k_3_6
    dydt[6*Lorg+5] = Jbb_core_6 * (A[-1]/V[-1]) + 1.0 * y[3*Lorg+2] * y[4*Lorg+3] * k_3_4- y[6*Lorg+5] * y[3*Lorg+2] * k_3_6

    #----component number 7, Trimer----
    dydt[6*Lorg+6] = Jbss_7 - Jssb_7 + 1.0 * y[2*Lorg+2] * y[5*Lorg+5] * k_3_6_surf
    dydt[6*Lorg+7] = (Jssb_7 - Jbss_7) * (A[0]/V[0]) + Jb1b2_7 * (A[1]/V[0]) + 1.0 * y[2*Lorg+3] * y[5*Lorg+6] * k_3_6
    for i in np.arange(1,Lorg-1):
        Jbb_minus1_7 = kbby_7[i] * (y[6*Lorg+7+i-1] - y[6*Lorg+7+i]) 
        Jbb_plus1_7 = kbby_7[i+1] * (y[6*Lorg+7+i+1] - y[6*Lorg+7+i]) 
        dydt[6*Lorg+7+i] = Jbb_minus1_7 * (A[i]/V[i]) + Jbb_plus1_7 * (A[i+1]/V[i]) + 1.0 * y[2*Lorg+3+i] * y[5*Lorg+6+i] * k_3_6
    dydt[7*Lorg+6] = Jbb_core_7 * (A[-1]/V[-1]) + 1.0 * y[3*Lorg+2] * y[6*Lorg+5] * k_3_6

    return dydt